The Dormand-Prince, Runge-Kutta integrator (5th order, with an embedded 4th order used for error estimation).
Appends the supplied solution point to the internal solution buffer.
Appends the supplied solution point to the internal solution buffer.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(inout) | :: | this |
The ode_integrator object. |
||
real(kind=real64), | intent(in) | :: | x |
The independent variable value. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | y |
The values of the dependent variables corresponding to x. |
|
class(errors), | intent(inout), | optional, | target | :: | err |
An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. Possible errors and warning messages that may be encountered are as follows.
|
Attempts an integration step for this integrator.
Attempts an integration step for this integrator.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(runge_kutta_45), | intent(inout) | :: | this |
The runge_kutta_45 object. |
||
class(ode_container), | intent(inout) | :: | sys |
The ode_container object containing the ODE's to integrate. |
||
real(kind=real64), | intent(in) | :: | h |
The current step size. |
||
real(kind=real64), | intent(in) | :: | x |
The current value of the independent variable. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | y |
An N-element array containing the current solution at x. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | f |
An N-element array containing the values of the derivatives at x. |
|
real(kind=real64), | intent(out), | dimension(:) | :: | yn |
An N-element array where this routine will write the next solution estimate at x + h. |
|
real(kind=real64), | intent(out), | dimension(:) | :: | fn |
An N-element array where this routine will write the next derivative estimate at x + h. |
|
real(kind=real64), | intent(out), | dimension(:) | :: | yerr |
An N-element array where this routine will write an estimate of the error in each equation. |
|
real(kind=real64), | intent(out), | dimension(:,:) | :: | k |
An N-by-NSTAGES matrix containing the derivatives at each stage. |
Clears the contents of the buffer.
Clears the contents of the buffer.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(inout) | :: | this |
The ode_integrator object. |
Computes the norm of the scaled error estimate.
Computes the norm of the scaled error estimate. A value less than one indicates a successful step. A value greater than one suggests that the results do not meet the requested tolerances.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(in) | :: | this |
The ode_integrator object. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | y |
The previously accepted solution array (N-element). |
|
real(kind=real64), | intent(in), | dimension(size(y)) | :: | yest |
An N-element array containing the next solution point estimate. |
|
real(kind=real64), | intent(in), | dimension(size(y)) | :: | yerr |
An N-element array containing the estimate of error for each equation. |
The norm of the scaled error.
Computes an estimate of an initial step size.
Computes an estimate of an initial step size.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(in) | :: | this |
The ode_integrator object. |
||
class(ode_container), | intent(in) | :: | sys |
The ode_container object containing the ODE's to integrate. |
||
real(kind=real64), | intent(in) | :: | xo |
The initial value of the independent variable. |
||
real(kind=real64), | intent(in) | :: | xf |
The final value of the independent variable. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | yo |
The initial values of the dependent variables (N-element). |
|
real(kind=real64), | intent(out), | dimension(size(yo)) | :: | fo |
An N-element array where the function values at xo will be written. |
|
real(kind=real64), | intent(out) | :: | h |
The initial step size estimate. |
Estimates the next step size.
Estimates the next step size based upon the current and previous error estimates.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(inout) | :: | this |
The ode_integrator object. |
||
real(kind=real64), | intent(in) | :: | e |
The norm of the current scaled error estimate. |
||
real(kind=real64), | intent(inout) | :: | eold |
The norm of the previous step's scaled error estimate. On output, this variable is updated. |
||
real(kind=real64), | intent(in) | :: | h |
The current step size. |
||
real(kind=real64), | intent(in) | :: | x |
The current independent variable value. |
||
class(errors), | intent(inout), | optional, | target | :: | err |
An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. Possible errors and warning messages that may be encountered are as follows.
|
The new step size estimate.
Gets the absolute error tolerance.
Gets the absolute error tolerance.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(in) | :: | this |
The ode_integrator object. |
The tolerance value.
Gets a value determining if the solver is allowed to overshoot the final value in the integration range.
Gets a value determining if the solver is allowed to overshoot the final value in the integration range.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(in) | :: | this |
The ode_integrator object. |
True if the solver can overshoot, and then interpolate to achieve the required final value; else, false thereby indicating the solver cannot overshoot.
Gets a logical parameter stating if this is a first-same-as-last (FSAL) integrator.
Gets a logical parameter stating if this is a first-same-as-last (FSAL) integrator.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(runge_kutta_45), | intent(in) | :: | this |
The runge_kutta_45 object. |
True for a FSAL integrator; else, false.
Gets the magnitude of the maximum allowed step size.
Gets the magnitude of the maximum allowed step size.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(in) | :: | this |
The ode_integrator object. |
The step size limit.
Gets the magnitude of the minimum allowed step size.
Gets the magnitude of the minimum allowed step size.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(in) | :: | this |
The ode_integrator object. |
The step size limit.
Gets the order of the integrator.
Gets the order of the integrator.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(runge_kutta_45), | intent(in) | :: | this |
The runge_kutta_45 object. |
The order.
Gets the relative error tolerance.
Gets the absolute error tolerance.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(in) | :: | this |
The ode_integrator object. |
The tolerance value.
Returns the solution computed by the integrator.
Returns the solution computed by the integrator stored as a matrix with the first column containing the values of the independent variable at which the solution was computed. The remaining columns contain the solutions for each of the integrated equations in the order in which they appear in the source routine. Notice, the solve routine must be called before this routine.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(in) | :: | this |
The ode_integrator object. |
The resulting solution matrix.
Gets the stage count for this integrator.
Gets the stage count for this integrator.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(runge_kutta_45), | intent(in) | :: | this |
The runge_kutta_45 object. |
The stage count.
Gets the limit on the number of integration steps.
Gets the limit on the number of integration steps that may be taken before the solver terminates.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(in) | :: | this |
The ode_integrator object. |
The step limit.
Gets the step size PI control parameter.
Gets the step size control parameter β used for PI control of the step size. A value of 0 provides a default step size controller (non-PI); however, a nonzero value of β provides PI control that improves stability, but comes with a potential for efficiency loss. A good estimate for a starting point for this parameter is β=0.4k where k is the order of the integrator.
The PI controller for step size is defined as follows. hn+1=fshne−αneβn−1
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(in) | :: | this |
The ode_integrator object. |
The control parameter.
Gets the step size safety factor.
Gets the safety factor (step size multiplier) used to provide a measure of control to the step size estimate such that hnew=fsh, where fs is this safety factor.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(in) | :: | this |
The ode_integrator object. |
The safety factor.
Performs the interpolation.
Performs the interpolation.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(runge_kutta_45), | intent(in) | :: | this |
The runge_kutta_45 object. |
||
real(kind=real64), | intent(in) | :: | x |
The value of the independent variable at which to compute the interpolation. |
||
real(kind=real64), | intent(in) | :: | xn |
The previous value of the independent variable at which the solution is computed. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | yn |
An N-element array containing the solution at xn. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | fn |
An N-element array containing the derivatives at xn. |
|
real(kind=real64), | intent(in) | :: | xn1 |
The value of the independent variable at xn + h. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | yn1 |
An N-element array containing the solution at xn + h. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | fn1 |
An N-element array containing the derivatives at xn + h. |
|
real(kind=real64), | intent(out), | dimension(:) | :: | y |
An N-element array where this routine will write the solution values interpolated at x. |
Sets up the interpolation process as the post-step action.
Sets up the interpolation process.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(runge_kutta_45), | intent(inout) | :: | this |
The runge_kutta_45 object. |
||
class(ode_container), | intent(inout) | :: | sys |
The ode_container object containing the ODE's to integrate. |
||
logical, | intent(in) | :: | dense |
Determines if dense output is requested (true); else, false. |
||
real(kind=real64), | intent(in) | :: | x |
The previous value of the independent variable. |
||
real(kind=real64), | intent(in) | :: | xn |
The current value of the independent variable. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | y |
An N-element array containing the solution at x. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | yn |
An N-element array containing the solution at xn. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | f |
An N-element array containing the derivatives at x. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | fn |
An N-element array containing the derivatives at xn. |
|
real(kind=real64), | intent(inout), | dimension(:,:) | :: | k |
An N-by-NSTAGES matrix containing the derivatives at each stage. |
Performs any pre-step actions.
Placeholder routine for any pre-step actions.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(runge_kutta_45), | intent(inout) | :: | this |
The runge_kutta_45 object. |
||
logical, | intent(in) | :: | prevs |
Defines the status of the previous step. The value is true if the previous step was successful; else, false if the previous step failed. |
||
class(ode_container), | intent(inout) | :: | sys |
The ode_container object containing the ODE's to integrate. |
||
real(kind=real64), | intent(in) | :: | h |
The current step size. |
||
real(kind=real64), | intent(in) | :: | x |
The current value of the independent variable. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | y |
An N-element array containing the current solution at x. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | f |
An N-element array containing the values of the derivatives at x. |
|
class(errors), | intent(inout), | optional, | target | :: | err |
An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. |
Sets the absolute error tolerance.
Sets the absolute error tolerance.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(inout) | :: | this |
The ode_integrator object. |
||
real(kind=real64), | intent(in) | :: | x |
The tolerance value. |
Sets a value determining if the solver is allowed to overshoot the final value in the integration range.
Sets a value determining if the solver is allowed to overshoot the final value in the integration range.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(inout) | :: | this |
The ode_integrator object. |
||
logical, | intent(in) | :: | x |
True if the solver can overshoot, and then interpolate to achieve the required final value; else, false thereby indicating the solver cannot overshoot. |
Sets the magnitude of the maximum allowed step size.
Sets the magnitude of the maximum allowed step size.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(inout) | :: | this |
The ode_integrator object. |
||
real(kind=real64), | intent(in) | :: | x |
The step size limit. |
Sets the magnitude of the minimum allowed step size.
Sets the magnitude of the minimum allowed step size.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(inout) | :: | this |
The ode_integrator object. |
||
real(kind=real64), | intent(in) | :: | x |
The step size limit. |
Sets the relative error tolerance.
Sets the absolute error tolerance.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(inout) | :: | this |
The ode_integrator object. |
||
real(kind=real64), | intent(in) | :: | x |
The tolerance value. |
Sets the limit on the number of integration steps.
Sets the limit on the number of integration steps that may be taken before the solver terminates.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(inout) | :: | this |
The ode_integrator object. |
||
integer(kind=int32), | intent(in) | :: | x |
The step limit. |
Sets the step size PI control parameter.
Sets the step size control parameter β used for PI control of the step size. A value of 0 provides a default step size controller (non-PI); however, a nonzero value of β provides PI control that improves stability, but comes with a potential for efficiency loss. A good estimate for a starting point for this parameter is β=0.4k where k is the order of the integrator.
The PI controller for step size is defined as follows. hn+1=fshne−αneβn−1
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(inout) | :: | this |
The ode_integrator object. |
||
real(kind=real64), | intent(in) | :: | x |
The control parameter |
Sets the step size safety factor.
Sets the safety factor (step size multiplier) used to provide a measure of control to the step size estimate such that hnew=fsh, where fs is this safety factor.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(inout) | :: | this |
The ode_integrator object. |
||
real(kind=real64), | intent(in) | :: | x |
The safety factor. |
Solves the supplied system of ODE's.
Solves the supplied system of ODE's.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(single_step_integrator), | intent(inout) | :: | this |
The single_step_integrator object. |
||
class(ode_container), | intent(inout) | :: | sys |
The ode_container object containing the ODE's to integrate. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | x |
An array of independent variable values at which to return the the solution to the ODE's. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | iv |
An array containing the initial values for each ODE. |
|
class(errors), | intent(inout), | optional, | target | :: | err |
An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. Possible errors and warning messages that may be encountered are as follows.
|